Top level modules

smuthi.simulation

Provide class to manage a simulation.

class smuthi.simulation.Simulation(layer_system=None, particle_list=None, initial_field=None, k_parallel='default', angular_resolution=0.008726646259971648, neff_waypoints=None, neff_imag=0.01, neff_max=None, neff_max_offset=1, neff_resolution=0.01, neff_minimal_branchpoint_distance=None, overwrite_default_contours=True, solver_type='LU', solver_tolerance=0.0001, store_coupling_matrix=True, coupling_matrix_lookup_resolution=None, coupling_matrix_interpolator_kind='linear', length_unit='length unit', input_file=None, output_dir='smuthi_output', save_after_run=False, log_to_file=False, log_to_terminal=True, check_circumscribing_spheres=True, do_sanity_check=True, periodicity=None, ewald_sum_separation_parameter='default', number_of_threads_periodic='default', use_pvwf_coupling=False, pvwf_coupling_neff_max=None, pvwf_coupling_neff_resolution=0.01)

Central class to manage a simulation.

Parameters:
  • layer_system (smuthi.layers.LayerSystem) – stratified medium
  • particle_list (list) – list of smuthi.particles.Particle objects
  • initial_field (smuthi.initial_field.InitialField) – initial field object
  • k_parallel (numpy.ndarray or str) – in-plane wavenumber for Sommerfeld integrals and field expansions. if ‘default’, use smuthi.fields.default_Sommerfeld_k_parallel_array
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
  • neff_waypoints (list or ndarray) – Used to set default k_parallel arrays. Corner points through which the contour runs This quantity is dimensionless (effective refractive index, will be multiplied by vacuum wavenumber) Multipole cut-off If not provided, reasonable waypoints are estimated.
  • neff_imag (float) – Used to set default k_parallel arrays. Extent of the contour into the negative imaginary direction (in terms of effective refractive index, n_eff=kappa/omega). Only needed when no neff_waypoints are provided
  • neff_max (float) – Used to set default k_parallel arrays. Truncation value of contour (in terms of effective refractive index). Only needed when no neff_waypoints are provided
  • neff_max_offset (float) – Used to set default k_parallel arrays. Use the last estimated singularity location plus this value (in terms of effective refractive index). Default=1 Only needed when no neff_waypoints are provided and if no value for neff_max is specified.
  • neff_resolution (float) – Used to set default k_parallel arrays. Resolution of contour, again in terms of effective refractive index
  • neff_minimal_branchpoint_distance (float) – Used to set default k_parallel arrays. Minimal distance that contour points shall have from branchpoint singularities (in terms of effective refractive index). This is only relevant if not deflected into imaginary. Default: One fifth of neff_resolution
  • overwrite_default_contours (bool) – If true (default), the default contours are written even if they have already been defined before
  • solver_type (str) – What solver type to use? Options: ‘LU’ for LU factorization, ‘gmres’ for GMRES iterative solver
  • coupling_matrix_lookup_resolution (float or None) – If type float, compute particle coupling by interpolation of a lookup table with that spacial resolution. If None (default), don’t use a lookup table but compute the coupling directly. This is more suitable for a small particle number.
  • coupling_matrix_interpolator_kind (str) – Set to ‘linear’ (default) or ‘cubic’ interpolation of the lookup table.
  • store_coupling_matrix (bool) – If True (default), the coupling matrix is stored. Otherwise it is recomputed on the fly during each iteration of the solver.
  • length_unit (str) – what is the physical length unit? has no influence on the computations
  • input_file (str) – path and filename of input file (for logging purposes)
  • output_dir (str) – path to folder where to export data
  • save_after_run (bool) – if true, the simulation object is exported to disc when over
  • log_to_file (bool) – if true, the simulation log will be written to a log file
  • log_to_terminal (bool) – if true, the simulation progress will be displayed in the terminal
  • check_circumscribing_spheres (bool) – if true, check all particles for overlapping circumscribing spheres and print a warning if detected
  • do_sanity_check (bool) – if true (default), check numerical input for some flaws. Warning: A passing sanity check does not guarantee correct numerical settings. For many particles, the sanity check might take some time and/or occupy large memory.
  • periodicity (tuple) – tuple (a1, a2) specifying two 3-dimensional lattice vectors in Carthesian coordinates with a1, a2 (numpy.ndarrays)
  • ewald_sum_separation_parameter (float) – Used to separate the real and reciprocal lattice sums to evaluate particle coupling in periodic lattices.
  • number_of_threads_periodic (int or str) – sets the number of threats used in a simulation with periodic particle arrangements if ‘default’, all available CPU cores are used if negative, all but number_of_threads_periodic are used
  • use_pvwf_coupling (bool) – If set to True, plane wave coupling is used to calculate the direct. Currently only possible in combination with direct solver strategy.
  • pvwf_coupling_neff_max (float) – Truncation neff for the integration contour of the PVWF coupling integral
  • pvwf_coupling_neff_resolution (float) – Discretization of the neff integral for PVWF coupling
circumscribing_spheres_disjoint()

Check if all circumscribing spheres are disjoint

initialize_linear_system()
largest_lateral_distance()

Compute the largest lateral distance between any two particles

print_simulation_header()
run()
sanity_check()

Check contour parameters for obvious problems

save(filename=None)

Export simulation object to disc.

Parameters:filename (str) – path and file name where to store data
set_default_Sommerfeld_contour()

Set the default Sommerfeld k_parallel array

set_default_angles()

Set the default polar and azimuthal angular arrays for pre-processing (i.e., initial field expansion)

set_default_contours()

Set the default initial field k_parallel array and the default Sommerfeld k_parallel array

set_default_initial_field_contour()

Set the default initial field k_parallel array

set_logging(log_to_terminal=None, log_to_file=None, log_filename=None)

Update logging behavior.

Parameters:
  • log_to_terminal (logical) – If true, print output to console.
  • log_to_file (logical) – If true, print output to file
  • log_filename (char) – If log_to_file is true, print output to a file with that name in the output directory. If the file already exists, it will be appended.

smuthi.initial_field

This module defines classes to represent the initial excitation.

class smuthi.initial_field.DipoleCollection(vacuum_wavelength, k_parallel_array='default', azimuthal_angles_array='default', angular_resolution=None, compute_swe_by_pwe=False, compute_dissipated_power_by_pwe=False)

Class for the representation of a set of point dipole sources. Use the append method to add DipoleSource objects.

Parameters:
  • vacuum_wavelength (float) – vacuum wavelength (length units)
  • k_parallel_array (numpy.ndarray or str) – In-plane wavenumber. If ‘default’, use smuthi.fields.default_initial_field_k_parallel_array
  • azimuthal_angles_array (numpy.ndarray or str) – Azimuthal angles for plane wave expansions If ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
  • compute_swe_by_pwe (bool) – If True, the initial field coefficients are computed through a plane wave expansion of the whole dipole collection field. This is slower for few dipoles and particles, but can become faster than the default for many dipoles and particles (default=False).
  • compute_dissipated_power_by_pwe (bool) – If True, evaluate dissipated power through a plane wave expansion of the whole scattered field. This is slower for few dipoles, but can be faster than the default for many dipoles (default=False).
append(dipole)

Add dipole to collection.

Parameters:dipole (DipoleSource) – Dipole object to add.
dissipated_power(particle_list, layer_system, k_parallel='default', azimuthal_angles='default', angular_resolution=None)

Compute the power that the dipole collection feeds into the system.

It is computed according to

\[P = \sum_i P_{0, i} + \frac{\omega}{2} \mathrm{Im} (\mathbf{\mu}_i^* \cdot \mathbf{E}_i(\mathbf{r}_i))\]

where \(P_{0,i}\) is the power that the i-th dipole would feed into an infinte homogeneous medium with the same refractive index as the layer that contains that dipole, \(\mathbf{r}_i\) is the location of the i-th dipole, \(\omega\) is the angular frequency, \(\mathbf{\mu}_i\) is the dipole moment and \(\mathbf{E}_i\) includes the reflections of the dipole field from the layer interfaces, as well as the scattered field from all particles and the fields from all other dipoles. In contrast to dissipated_power_alternative, this routine uses the particle coupling routines and might be faster for many particles and few dipoles.

Parameters:
  • particle_list (list of smuthi.particles.Particle objects) – scattering particles
  • layer_system (smuthi.layers.LayerSystem) – stratified medium
  • k_parallel (ndarray or str) – array of in-plane wavenumbers for plane wave expansions. If ‘default’, use smuthi.fields.default_initial_field_k_parallel_array
  • azimuthal_angles (ndarray or str) – array of azimuthal angles for plane wave expansions. If ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
Returns:

dissipated power of each dipole (list of floats)

dissipated_power_alternative(particle_list, layer_system, k_parallel='default', azimuthal_angles='default', angular_resolution=None)

Compute the power that the dipole collection feeds into the system.

It is computed according to

\[P = \sum_i P_{0, i} + \frac{\omega}{2} \mathrm{Im} (\mathbf{\mu}_i^* \cdot \mathbf{E}_i(\mathbf{r}_i))\]

where \(P_{0,i}\) is the power that the i-th dipole would feed into an infinte homogeneous medium with the same refractive index as the layer that contains that dipole, \(\mathbf{r}_i\) is the location of the i-th dipole, \(\omega\) is the angular frequency, \(\mathbf{\mu}_i\) is the dipole moment and \(\mathbf{E}_i\) includes the reflections of the dipole field from the layer interfaces, as well as the scattered field from all particles and the fields from all other dipoles. In contrast to dissipated_power, this routine uses the scattered field piecewise expansion and might be faster for few particles or many dipoles.

Parameters:
  • particle_list (list of smuthi.particles.Particle objects) – scattering particles
  • layer_system (smuthi.layers.LayerSystem) – stratified medium
  • k_parallel (ndarray or str) – array of in-plane wavenumbers for plane wave expansions. If ‘default’, use smuthi.fields.default_initial_field_k_parallel_array
  • azimuthal_angles (ndarray or str) – array of azimuthal angles for plane wave expansions. If ‘default’, use smuthi.fields.default_azimuthal_angles
  • angular_resolution (float) – If provided, angular arrays are generated with this angular resolution over the default angular range
Returns:

dissipated power of each dipole (list of floats)

electric_field(x, y, z, layer_system)

Evaluate the complex electric field of the dipole collection.

Parameters:
  • x (array like) – Array of x-values where to evaluate the field (length unit)
  • y (array like) – Array of y-values where to evaluate the field (length unit)
  • z (array like) – Array of z-values where to evaluate the field (length unit)
  • layer_system (smuthi.layer.LayerSystem) – Stratified medium
Returns
Tuple (E_x, E_y, E_z) of electric field values
magnetic_field(x, y, z, layer_system)

Evaluate the complex magnetic field of the dipole collection.

Parameters:
  • x (array like) – Array of x-values where to evaluate the field (length unit)
  • y (array like) – Array of y-values where to evaluate the field (length unit)
  • z (array like) – Array of z-values where to evaluate the field (length unit)
  • layer_system (smuthi.layer.LayerSystem) – Stratified medium
Returns
Tuple (H_x, H_y, H_z) of magnetic field values
piecewise_field_expansion(layer_system)

Compute a piecewise field expansion of the dipole collection..

Parameters:layer_system (smuthi.layer.LayerSystem) – stratified medium
Returns:smuthi.field_expansion.PiecewiseWaveExpansion object
plane_wave_expansion = functools.partial(<bound method Memoize.__call__ of <smuthi.utility.memoizing.Memoize object>>, None)
spherical_wave_expansion(particle, layer_system)

Regular spherical wave expansion of the dipole collection including layer system response, at the locations of the particles. If self.compute_swe_by_pwe is True, use the dipole collection plane wave expansion, otherwise use the individual dipoles spherical_wave_expansion method.

Parameters:
  • particle (smuthi.particles.Particle) – particle relative to which the swe is computed
  • layer_system (smuthi.layer.LayerSystem) – stratified medium
Returns:

regular smuthi.field_expansion.SphericalWaveExpansion object

class smuthi.initial_field.DipoleSource(vacuum_wavelength, dipole_moment, position, k_parallel_array='default', azimuthal_angles_array='default')

Class for the representation of a single point dipole source.

Parameters:
  • vacuum_wavelength (float) – vacuum wavelength (length units)
  • dipole_moment (list or tuple) – (x, y, z)-coordinates of dipole moment vector
  • position (list or tuple) – (x, y, z)-coordinates of dipole position
  • k_parallel_array (numpy.ndarray or str) – In-plane wavenumber. If ‘default’, use smuthi.fields.default_initial_field_k_parallel_array
  • azimuthal_angles_array (numpy.ndarray or str) – Azimuthal angles for plane wave expansions If ‘default’, use smuthi.fields.default_azimuthal_angles
check_dissipated_power_homogeneous_background(layer_system)
current()

The current density takes the form

\[\mathbf{j}(\mathbf{r}) = \delta(\mathbf{r} - \mathbf{r}_D) \mathbf{j}_D,\]

where \(\mathbf{j}_D = -j \omega \mathbf{\mu}\), \(\mathbf{r}_D\) is the location of the dipole, \(\omega\) is the angular frequency and \(\mathbf{\mu}\) is the dipole moment. For further details, see ‘Principles of nano optics’ by Novotny and Hecht.

Returns:List of [x, y, z]-components of current density vector \(\mathbf{j}_D\)
dissipated_power(particle_list, layer_system, show_progress=True)

Compute the power that the dipole feeds into the system.

It is computed according to

\[P = P_0 + \frac{\omega}{2} \mathrm{Im} (\mathbf{\mu}^* \cdot \mathbf{E}(\mathbf{r}_D))\]

where \(P_0\) is the power that the dipole would feed into an infinte homogeneous medium with the same refractive index as the layer that contains the dipole, \(\mathbf{r}_D\) is the location of the dipole, \(\omega\) is the angular frequency, \(\mathbf{\mu}\) is the dipole moment and \(\mathbf{E}\) includes the reflections of the dipole field from the layer interfaces, as well as the scattered field from all particles.

Parameters:
  • particle_list (list of smuthi.particles.Particle objects) – scattering particles
  • layer_system (smuthi.layers.LayerSystem) – stratified medium
  • show_progress (bool) – if true, display progress
Returns:

dissipated power as float

dissipated_power_alternative(particle_list, layer_system)

Compute the power that the dipole feeds into the system.

It is computed according to

\[P = P_0 + \frac{\omega}{2} \mathrm{Im} (\mathbf{\mu}^* \cdot \mathbf{E}(\mathbf{r}_D))\]

where \(P_0\) is the power that the dipole would feed into an infinte homogeneous medium with the same refractive index as the layer that contains the dipole, \(\mathbf{r}_D\) is the location of the dipole, \(\omega\) is the angular frequency, \(\mathbf{\mu}\) is the dipole moment and \(\mathbf{E}\) includes the reflections of the dipole field from the layer interfaces, as well as the scattered field from all particles. In contrast to dissipated_power, this routine relies on the scattered field piecewise expansion and and might thus be slower.

Parameters:
  • particle_list (list of smuthi.particles.Particle objects) – scattering particles
  • layer_system (smuthi.layers.LayerSystem) – stratified medium
Returns:

dissipated power as float

dissipated_power_homogeneous_background(layer_system)

Compute the power that the dipole would radiate in an infinite homogeneous medium of the same refractive index as the layer that contains the dipole.

\[P_0 = \frac{|\mathbf{\mu}| k \omega^3}{12 \pi}\]
Parameters:layer_system (smuthi.layers.LayerSystem) – stratified medium
Returns:power (float)
electric_field(x, y, z, layer_system, include_direct_field=True, include_layer_response=True)

Evaluate the complex electric field of the dipole source.

Parameters:
  • x (array like) – Array of x-values where to evaluate the field (length unit)
  • y (array like) – Array of y-values where to evaluate the field (length unit)
  • z (array like) – Array of z-values where to evaluate the field (length unit)
  • layer_system (smuthi.layer.LayerSystem) – Stratified medium
  • include_direct_field (bool) – if True (default), the direct dipole field is included. otherwise, only the layer response of the dipole field is returned.
  • include_layer_response (bool) – if True (default), the layer response of the dipole field is included. otherwise, only the direct dipole field is returned.
Returns
Tuple (E_x, E_y, E_z) of electric field values
magnetic_field(x, y, z, layer_system, include_direct_field=True, include_layer_response=True)

Evaluate the complex magnetic field of the dipole source.

Parameters:
  • x (array like) – Array of x-values where to evaluate the field (length unit)
  • y (array like) – Array of y-values where to evaluate the field (length unit)
  • z (array like) – Array of z-values where to evaluate the field (length unit)
  • layer_system (smuthi.layer.LayerSystem) – Stratified medium
  • include_direct_field (bool) – if True (default), the direct dipole field is included. otherwise, only the layer response of the dipole field is returned.
  • include_layer_response (bool) – if True (default), the layer response of the dipole field is included. otherwise, only the direct dipole field is returned.
Returns
Tuple (H_x, H_y, H_z) of electric field values
outgoing_spherical_wave_expansion(layer_system)

The dipole field as an expansion in spherical vector wave functions.

Parameters:layer_system (smuthi.layers.LayerSystem) – stratified medium
Returns:outgoing smuthi.field_expansion.SphericalWaveExpansion object
piecewise_field_expansion(layer_system, include_direct_field=True, include_layer_response=True)

Compute a piecewise field expansion of the dipole field.

Parameters:
  • layer_system (smuthi.layer.LayerSystem) – stratified medium
  • include_direct_field (bool) – if True (default), the direct dipole field is included. otherwise, only the layer response of the dipole field is returned.
  • include_layer_response (bool) – if True (default), the layer response of the dipole field is included. otherwise, only the direct dipole field is returned.
Returns:

smuthi.field_expansion.PiecewiseWaveExpansion object

plane_wave_expansion(layer_system, i, k_parallel_array=None, azimuthal_angles_array=None)

Plane wave expansion of the dipole field.

Parameters:
  • layer_system (smuthi.layer.LayerSystem) – stratified medium
  • i (int) – layer number in which to evaluate the expansion
  • k_parallel_array (numpy.ndarray) – in-plane wavenumber array for the expansion. if none specified, self.k_parallel_array is used
  • azimuthal_angles_array (numpy.ndarray) – azimuthal angles for the expansion. if none specified, self.azimuthal_angles_array is used
Returns:

tuple of to smuthi.field_expansion.PlaneWaveExpansion objects, one for upgoing and one for downgoing component

spherical_wave_expansion(particle, layer_system)

Regular spherical wave expansion of the wave including layer system response, at the locations of the particles.

Parameters:
  • particle (smuthi.particles.Particle) – particle relative to which the swe is computed
  • layer_system (smuthi.layer.LayerSystem) – stratified medium
Returns:

regular smuthi.field_expansion.SphericalWaveExpansion object

class smuthi.initial_field.GaussianBeam(vacuum_wavelength, polar_angle, azimuthal_angle, polarization, beam_waist, k_parallel_array='default', azimuthal_angles_array='default', amplitude=1, reference_point=None)

Class for the representation of a Gaussian beam as initial field.

initial_intensity(layer_system)

Evaluate the incoming intensity of the initial field.

Parameters:layer_system (smuthi.layers.LayerSystem) – Stratified medium
Returns:A smuthi.field_expansion.FarField object holding the initial intensity information.
plane_wave_expansion(layer_system, i, k_parallel_array=None, azimuthal_angles_array=None)

Plane wave expansion of the Gaussian beam.

Parameters:
  • layer_system (smuthi.layer.LayerSystem) – stratified medium
  • i (int) – layer number in which to evaluate the expansion
  • k_parallel_array (numpy.ndarray) – in-plane wavenumber array for the expansion. if none specified, self.k_parallel_array is used
  • azimuthal_angles_array (numpy.ndarray) – azimuthal angles for the expansion. if none specified, self.azimuthal_angles_array is used
Returns:

tuple of to smuthi.field_expansion.PlaneWaveExpansion objects, one for upgoing and one for downgoing component

propagated_far_field(layer_system)

Evaluate the far field intensity of the reflected / transmitted initial field.

Parameters:layer_system (smuthi.layers.LayerSystem) – Stratified medium
Returns:A tuple of smuthi.field_expansion.FarField objects, one for forward (i.e., into the top hemisphere) and one for backward propagation (bottom hemisphere).
class smuthi.initial_field.InitialField(vacuum_wavelength)

Base class for initial field classes

angular_frequency()

Angular frequency.

Returns:Angular frequency (float) according to the vacuum wavelength in units of c=1.
get_azimuthal_angles_array()

Get azimuthal angles array which is either the default array or the one stored in the object

get_k_parallel_array()

Get k_parallel array which is either the default array or the one stored in the object

piecewise_field_expansion(layer_system)
plane_wave_expansion(layer_system, i)

Virtual method to be overwritten.

spherical_wave_expansion(particle, layer_system)

Virtual method to be overwritten.

class smuthi.initial_field.InitialPropagatingWave(vacuum_wavelength, polar_angle, azimuthal_angle, polarization, amplitude=1, reference_point=None)

Base class for plane waves and Gaussian beams

Parameters:
  • vacuum_wavelength (float) –
  • polar_angle (float) – polar propagation angle (0 means, parallel to z-axis)
  • azimuthal_angle (float) – azimuthal propagation angle (0 means, in x-z plane)
  • polarization (int) – 0 for TE/s, 1 for TM/p
  • amplitude (float or complex) – Electric field amplitude
  • reference_point (list) – Location where electric field of incoming wave equals amplitude
electric_field(x, y, z, layer_system)

Evaluate the complex electric field corresponding to the wave.

Parameters:
  • x (array like) – Array of x-values where to evaluate the field (length unit)
  • y (array like) – Array of y-values where to evaluate the field (length unit)
  • z (array like) – Array of z-values where to evaluate the field (length unit)
  • layer_system (smuthi.layer.LayerSystem) – Stratified medium
Returns
Tuple (E_x, E_y, E_z) of electric field values
magnetic_field(x, y, z, layer_system)

Evaluate the complex magnetic field corresponding to the wave.

Parameters:
  • x (array like) – Array of x-values where to evaluate the field (length unit)
  • y (array like) – Array of y-values where to evaluate the field (length unit)
  • z (array like) – Array of z-values where to evaluate the field (length unit)
  • layer_system (smuthi.layer.LayerSystem) – Stratified medium
Returns
Tuple (H_x, H_y, H_z) of magnetic field values
piecewise_field_expansion(layer_system)

Compute a piecewise field expansion of the initial field.

Parameters:layer_system (smuthi.layer.LayerSystem) – stratified medium
Returns:smuthi.field_expansion.PiecewiseWaveExpansion object
spherical_wave_expansion(particle, layer_system)

Regular spherical wave expansion of the wave including layer system response, at the locations of the particles.

Parameters:
  • particle (smuthi.particles.Particle) – particle relative to which the swe is computed
  • layer_system (smuthi.layer.LayerSystem) – stratified medium
Returns:

regular smuthi.field_expansion.SphericalWaveExpansion object

class smuthi.initial_field.PlaneWave(vacuum_wavelength, polar_angle, azimuthal_angle, polarization, amplitude=1, reference_point=None)

Class for the representation of a plane wave as initial field.

Parameters:
  • vacuum_wavelength (float) –
  • polar_angle (float) – polar angle of k-vector (0 means, k is parallel to z-axis)
  • azimuthal_angle (float) – azimuthal angle of k-vector (0 means, k is in x-z plane)
  • polarization (int) – 0 for TE/s, 1 for TM/p
  • amplitude (float or complex) – Plane wave amplitude at reference point
  • reference_point (list) – Location where electric field of incoming wave equals amplitude
plane_wave_expansion(layer_system, i)

Plane wave expansion for the plane wave including its layer system response. As it already is a plane wave, the plane wave expansion is somehow trivial (containing only one partial wave, i.e., a discrete plane wave expansion).

Parameters:
  • layer_system (smuthi.layers.LayerSystem) – Layer system object
  • i (int) – layer number in which the plane wave expansion is valid
Returns:

Tuple of smuthi.field_expansion.PlaneWaveExpansion objects. The first element is an upgoing PWE, whereas the second element is a downgoing PWE.

smuthi.layers

Provide class for the representation of planar layer systems.

class smuthi.layers.LayerSystem(thicknesses=None, refractive_indices=None)

Stack of planar layers.

Parameters:
  • thicknesses (list) – layer thicknesses, first and last are semi inf and set to 0 (length unit)
  • refractive_indices (list) – complex refractive indices in the form n+jk
is_degenerate()

Returns True if the layer system consists of only two layers of the same material. This function is useful for detecting if layer mediated coupling can be omitted when calculating the coupling between particles.

Returns:True if layer system is degenerate. False otherwise.
layer_number(z)

Return number of layer that contains point [0,0,z]

If z is on the interface, the higher layer number is selected.

Parameters:z (float) – z-coordinate of query point (length unit)
Returns:number of layer containing z
lower_zlimit(i)

Return the z-coordinate of lower boundary

The coordinate system is defined such that z=0 corresponds to the interface between layer 0 and layer 1.

Parameters:i (int) – index of layer in question (must be between 0 and number_of_layers-1)
Returns:z-coordinate of lower boundary
number_of_layers()

Return total number of layers

Returns:number of layers
reference_z(i)

Return the anchor point’s z-coordinate.

The coordinate system is defined such that z=0 corresponds to the interface between layer 0 and layer 1.

Parameters:i (int) – index of layer in question (must be between 0 and number_of_layers-1)
Returns:anchor point’s z-coordinate
response(pwe, from_layer, to_layer)

Evaluate the layer system response to an electromagnetic excitation inside the layer system.

Parameters:
  • pwe (tuple or smuthi.field_expansion.PlaneWaveExpansion) – Either specify a PlaneWaveExpansion object that that represents the electromagnetic excitation, or a tuple of two PlaneWaveExpansion objects representing the upwards- and downwards propagating partial waves of the excitation.
  • from_layer (int) – Layer number in which the excitation is located
  • to_layer (int) – Layer number in which the layer response is to be evaluated
Returns:

Tuple (pwe_up, pwe_sown) of PlaneWaveExpansion objects representing the layer system response to the excitation.

upper_zlimit(i)

Return the z-coordinate of upper boundary.

The coordinate system is defined such that z=0 corresponds to the interface between layer 0 and layer 1.

Parameters:i (int) – index of layer in question (must be between 0 and number_of_layers-1)
Returns:z-coordinate of upper boundary
wavenumber(layer_number, vacuum_wavelength)
Parameters:
  • layer_number (int) – number of layer in question
  • vacuum_wavelength (float) – vacuum wavelength
Returns:

wavenumber in that layer as float

smuthi.layers.fresnel_r(pol, kz1, kz2, n1, n2)

Fresnel reflection coefficient.

Parameters:
  • pol (int) – polarization (0=TE, 1=TM)
  • kz1 (float or array) – incoming wave’s z-wavenumber (k*cos(alpha1))
  • kz2 (float or array) – transmitted wave’s z-wavenumber (k*cos(alpha2))
  • n1 (float or complex) – first medium’s complex refractive index (n+ik)
  • n2 (float or complex) – second medium’s complex refractive index (n+ik)
Returns:

Complex Fresnel reflection coefficient (float or array)

smuthi.layers.fresnel_t(pol, kz1, kz2, n1, n2)

Fresnel transmission coefficient.

Parameters:
  • pol (int) – polarization (0=TE, 1=TM)
  • kz1 (float or array) – incoming wave’s z-wavenumber (k*cos(alpha1))
  • kz2 (float or array) – transmitted wave’s z-wavenumber (k*cos(alpha2))
  • n1 (float or complex) – first medium’s complex refractive index (n+ik)
  • n2 (float or complex) – second medium’s complex refractive index (n+ik)
Returns:

Complex Fresnel transmission coefficient (float or array)

smuthi.layers.interface_transition_matrix(pol, kz1, kz2, n1, n2)

Interface transition matrix to be used in the Transfer matrix algorithm.

Parameters:
  • pol (int) – polarization (0=TE, 1=TM)
  • kz1 (float or array) – incoming wave’s z-wavenumber (k*cos(alpha1))
  • kz2 (float or array) – transmitted wave’s z-wavenumber (k*cos(alpha2))
  • n1 (float or complex) – first medium’s complex refractive index (n+ik)
  • n2 (float or complex) – second medium’s complex refractive index (n+ik)
Returns:

Interface transition matrix as 2x2 numpy array or as 2x2 mpmath.matrix

smuthi.layers.layer_propagation_matrix(kz, d)

Layer propagation matrix to be used in the Transfer matrix algorithm.

Parameters:
  • kz (float or complex) – z-wavenumber (k*cos(alpha))
  • d (float) – thickness of layer
Returns:

Layer propagation matrix as 2x2 numpy array or as 2x2 mpmath.matrix

smuthi.layers.layersystem_scattering_matrix(pol, layer_d, layer_n, kpar, omega)

Scattering matrix of a planarly layered medium.

Parameters:
  • pol (int) – polarization(0=TE, 1=TM)
  • layer_d (list) – layer thicknesses
  • layer_n (list) – complex layer refractive indices
  • kpar (float) – in-plane wavenumber
  • omega (float) – angular frequency in units of c=1: omega=2*pi/lambda
Returns:

Scattering matrix as 2x2 numpy array or as 2x2 mpmath.matrix

smuthi.layers.layersystem_transfer_matrix(pol, layer_d, layer_n, kpar, omega)

Transfer matrix of a planarly layered medium.

Parameters:
  • pol (int) – polarization(0=TE, 1=TM)
  • layer_d (list) – layer thicknesses
  • layer_n (list) – complex layer refractive indices
  • kpar (float) – in-plane wavenumber
  • omega (float) – angular frequency in units of c=1: omega=2*pi/lambda
Returns:

Transfer matrix as 2x2 numpy array or as 2x2 mpmath.matrix

smuthi.layers.matrix_inverse(m)
Parameters:m (mpmath.matrix or numpy.ndarray) – matrix to invert
Returns:inverse of m with same data type as m1 and m2
smuthi.layers.matrix_product(m1, m2)
Parameters:
  • m1 (mpmath.matrix or numpy.ndarray) – first matrix
  • m2 (mpmath.matrix or numpy.ndarray) – second matrix
Returns:

matrix product m1 * m2 with same data type as m1 and m2

smuthi.layers.set_precision(prec=None)

Set the numerical precision of the layer system response. You can use this to evaluate the layer response of unstable systems, for example in the case of evanescent waves in very thick layers. Calculations take longer time if the precision is set to a value other than None (default).

Parameters:prec (None or int) – If None, calculations are done using standard double precision. If int, that many decimal digits are considered in the calculations, using the mpmath package.

smuthi.particles

Classes for the representation of scattering particles.

class smuthi.particles.AnisotropicSphere(position=None, euler_angles=None, polar_angle=0, azimuthal_angle=0, refractive_index=(1+0j), radius=1, refractive_index_z=(2+0j), l_max=None, m_max=None, n_rank=None)

Particle subclass for anisotropic spheres.

Parameters:
  • position (list) – Particle position in the format [x, y, z] (length unit)
  • euler_angles (list) – Euler angles [alpha, beta, gamma] in (zy’z’’-convention) in radian. Alternatively, you can specify the polar and azimuthal angle of the axis of revolution.
  • polar_angle (float) – Polar angle of axis of revolution.
  • azimuthal_angle (float) – Azimuthal angle of axis of revolution.
  • refractive_index (complex) – Complex refractive index of particle in x-y plane (if not rotated)
  • refractive_index_z (complex) – Complex refractive index of particle along z-axis (if not rotated)
  • radius (float) – Sphere radius
  • l_max (int) – Maximal multipole degree used for the spherical wave expansion of incoming and scattered field
  • m_max (int) – Maximal multipole order used for the spherical wave expansion of incoming and scattered field
  • n_rank (int) – Maximal multipole order used for in NFMDS (default: l_max + 5)
circumscribing_sphere_radius()

Virtual method to be overwritten

compute_t_matrix(vacuum_wavelength, n_medium)

Return the T-matrix of a particle.

Parameters:
  • vacuum_wavelength (float) –
  • n_medium (float or complex) – Refractive index of surrounding medium
  • particle (smuthi.particles.Particle) – Particle object
Returns:

T-matrix as ndarray

class smuthi.particles.CustomParticle(position=None, euler_angles=None, polar_angle=0, azimuthal_angle=0, refractive_index=(1+0j), geometry_filename=None, scale=1, l_max=None, m_max=None, n_rank=None)

Particle subclass for custom particle shapes defined via FEM file.

Parameters:
  • position (list) – Particle position in the format [x, y, z] (length unit)
  • euler_angles (list) – Euler angles [alpha, beta, gamma] in (zy’z’’-convention) in radian. Alternatively, you can specify the polar and azimuthal angle of the axis of revolution.
  • polar_angle (float) – Polar angle of axis of revolution.
  • azimuthal_angle (float) – Azimuthal angle of axis of revolution.
  • geometry_filename (string) – Path to FEM file
  • scale (float) – Scaling factor for particle dimensions (relative to provided dimensions)
  • l_max (int) – Maximal multipole degree used for the spherical wave expansion of incoming and scattered field
  • m_max (int) – Maximal multipole order used for the spherical wave expansion of incoming and scattered field
  • n_rank (int) – Maximal multipole order used for in NFMDS (default: l_max + 5)
circumscribing_sphere_radius()

Virtual method to be overwritten

compute_t_matrix(vacuum_wavelength, n_medium)

Return the T-matrix of a particle.

Parameters:
  • vacuum_wavelength (float) –
  • n_medium (float or complex) – Refractive index of surrounding medium
  • particle (smuthi.particles.Particle) – Particle object
Returns:

T-matrix as ndarray

class smuthi.particles.FiniteCylinder(position=None, euler_angles=None, polar_angle=0, azimuthal_angle=0, refractive_index=(1+0j), cylinder_radius=1, cylinder_height=1, l_max=None, m_max=None, n_rank=None, use_python_tmatrix=False, nint=100)

Particle subclass for finite cylinders.

Parameters:
  • position (list) – Particle position in the format [x, y, z] (length unit)
  • euler_angles (list) – Euler angles [alpha, beta, gamma] in (zy’z’’-convention) in radian. Alternatively, you can specify the polar and azimuthal angle of the axis of revolution.
  • polar_angle (float) – Polar angle of axis of revolution.
  • azimuthal_angle (float) – Azimuthal angle of axis of revolution.
  • refractive_index (complex) – Complex refractive index of particle
  • cylinder_radius (float) – Radius of cylinder (length unit)
  • cylinder_height (float) – Height of cylinder, in z-direction if not rotated (length unit)
  • l_max (int) – Maximal multipole degree used for the spherical wave expansion of incoming and scattered field
  • m_max (int) – Maximal multipole order used for the spherical wave expansion of incoming and scattered field
  • n_rank (int) – Maximal multipole order used for in NFMDS (default: l_max + 5)
  • use_python_tmatrix (bool) – If true, use Alan Zhan’s Python code to compute the T-matrix rather than NFM-DS
  • nint (int) – Number of angles used in integral (only for python t-mnatrix)
circumscribing_sphere_radius()

Virtual method to be overwritten

compute_t_matrix(vacuum_wavelength, n_medium)

Return the T-matrix of a particle.

Parameters:
  • vacuum_wavelength (float) –
  • n_medium (float or complex) – Refractive index of surrounding medium
  • particle (smuthi.particles.Particle) – Particle object
Returns:

T-matrix as ndarray

class smuthi.particles.LayeredSpheroid(position=None, euler_angles=None, polar_angle=0, azimuthal_angle=0, layer_refractive_indices=(1+0j), layer_semi_axes_c=1, layer_semi_axes_a=1, l_max=None, m_max=None, n_rank=None)

Particle subclass for layered spheroid.

Parameters:
  • position (list) – Particle position in the format [x, y, z] (length unit)
  • euler_angles (list) – Euler angles [alpha, beta, gamma] in (zy’z’’-convention) in radian. Alternatively, you can specify the polar and azimuthal angle of the axis of revolution.
  • polar_angle (float) – Polar angle of axis of revolution.
  • azimuthal_angle (float) – Azimuthal angle of axis of revolution.
  • layer_refractive_indices (complex) – Complex refractive index of particle
  • layer_semi_axes_c (float) – Spheroid half axis in direction of axis of revolution (z-axis if not rotated)
  • layer_semi_axes_a (float) – Spheroid half axis in lateral direction (x- and y-axis if not rotated)
  • l_max (int) – Maximal multipole degree used for the spherical wave expansion of incoming and scattered field
  • m_max (int) – Maximal multipole order used for the spherical wave expansion of incoming and scattered field
  • n_rank (int) – Maximal multipole order used in NFMDS (default: l_max + 5)
compute_t_matrix(vacuum_wavelength, n_medium)

Return the T-matrix of a particle.

Parameters:
  • vacuum_wavelength (float) –
  • n_medium (float or complex) – Refractive index of surrounding medium
  • particle (smuthi.particles.Particle) – Particle object
Returns:

T-matrix as ndarray

class smuthi.particles.Particle(position=None, euler_angles=None, refractive_index=(1+0j), l_max=None, m_max=None)

Base class for scattering particles.

Parameters:
  • position (list) – Particle position in the format [x, y, z] (length unit)
  • euler_angles (list) – Particle Euler angles in the format [alpha, beta, gamma]
  • refractive_index (complex) – Complex refractive index of particle
  • l_max (int) – Maximal multipole degree used for the spherical wave expansion of incoming and scattered field
  • m_max (int) – Maximal multipole order used for the spherical wave expansion of incoming and scattered field
circumscribing_sphere_radius()

Virtual method to be overwritten

compute_t_matrix(vacuum_wavelength, n_medium)

Return the T-matrix of a particle.

Parameters:
  • vacuum_wavelength (float) –
  • n_medium (float or complex) – Refractive index of surrounding medium
  • particle (smuthi.particles.Particle) – Particle object
Returns:

T-matrix as ndarray

is_inside(x, y, z)

Virtual method to be overwritten. Until all child classes implement it: return False

is_outside(x, y, z)

Virtual method to be overwritten. Until all child classes implement it: return True

class smuthi.particles.Sphere(position=None, refractive_index=(1+0j), radius=1, l_max=None, m_max=None)

Particle subclass for spheres.

Parameters:
  • position (list) – Particle position in the format [x, y, z] (length unit)
  • refractive_index (complex) – Complex refractive index of particle
  • radius (float) – Particle radius (length unit)
  • l_max (int) – Maximal multipole degree used for the spherical wave expansion of incoming and scattered field
  • m_max (int) – Maximal multipole order used for the spherical wave expansion of incoming and scattered field
circumscribing_sphere_radius()

Virtual method to be overwritten

compute_t_matrix(vacuum_wavelength, n_medium)

Return the T-matrix of a particle.

Parameters:
  • vacuum_wavelength (float) –
  • n_medium (float or complex) – Refractive index of surrounding medium
  • particle (smuthi.particles.Particle) – Particle object
Returns:

T-matrix as ndarray

is_inside(x, y, z)

Virtual method to be overwritten. Until all child classes implement it: return False

is_outside(x, y, z)

Virtual method to be overwritten. Until all child classes implement it: return True

class smuthi.particles.Spheroid(position=None, euler_angles=None, polar_angle=0, azimuthal_angle=0, refractive_index=(1+0j), semi_axis_c=1, semi_axis_a=1, l_max=None, m_max=None, n_rank=None)

Particle subclass for spheroids.

Parameters:
  • position (list) – Particle position in the format [x, y, z] (length unit)
  • euler_angles (list) – Euler angles [alpha, beta, gamma] in (zy’z’’-convention) in radian. Alternatively, you can specify the polar and azimuthal angle of the axis of revolution.
  • polar_angle (float) – Polar angle of axis of revolution.
  • azimuthal_angle (float) – Azimuthal angle of axis of revolution.
  • refractive_index (complex) – Complex refractive index of particle
  • semi_axis_c (float) – Spheroid half axis in direction of axis of revolution (z-axis if not rotated)
  • semi_axis_a (float) – Spheroid half axis in lateral direction (x- and y-axis if not rotated)
  • l_max (int) – Maximal multipole degree used for the spherical wave expansion of incoming and scattered field
  • m_max (int) – Maximal multipole order used for the spherical wave expansion of incoming and scattered field
  • n_rank (int) – Maximal multipole order used in NFMDS (default: l_max + 5)
circumscribing_sphere_radius()

Virtual method to be overwritten

compute_t_matrix(vacuum_wavelength, n_medium)

Return the T-matrix of a particle.

Parameters:
  • vacuum_wavelength (float) –
  • n_medium (float or complex) – Refractive index of surrounding medium
  • particle (smuthi.particles.Particle) – Particle object
Returns:

T-matrix as ndarray